{
 "cells": [
  {
   "cell_type": "markdown",
   "metadata": {},
   "source": [
    "## Using Variational Estimates to Initialize the NUTS-HMC Sampler\n",
    "\n",
    "In this example we show how to use the parameter estimates return by Stan's variational inference algorithm\n",
    "as the initial parameter values for Stan's NUTS-HMC sampler.\n",
    "By default, the sampler algorithm randomly initializes all model parameters in the range uniform\\[-2, 2\\].  When the true parameter value is outside of this range, starting from the ADVI estimates will speed up and improve adaptation.\n",
    "\n",
    "### Model and data\n",
    "\n",
    "The Stan model and data are taken from the [posteriordb package](https://github.com/stan-dev/posteriordb).\n",
    "\n",
    "We use the [blr model](https://github.com/stan-dev/posteriordb/blob/master/posterior_database/models/stan/blr.stan),\n",
    "a Bayesian standard linear regression model with noninformative priors,\n",
    "and its corresponding simulated dataset [sblri.json](https://github.com/stan-dev/posteriordb/blob/master/posterior_database/data/data/sblri.json.zip),\n",
    "which was simulated via script [sblr.R](https://github.com/stan-dev/posteriordb/blob/master/posterior_database/data/data-raw/sblr/sblr.R).\n",
    "For conveince, we have copied the posteriordb model and data to this directory, in files `blr.stan` and `sblri.json`."
   ]
  },
  {
   "cell_type": "code",
   "execution_count": null,
   "metadata": {},
   "outputs": [],
   "source": [
    "import os\n",
    "from cmdstanpy import CmdStanModel\n",
    "\n",
    "stan_file = 'blr.stan' # basic linear regression\n",
    "data_file = 'sblri.json' # simulated data\n",
    "\n",
    "model = CmdStanModel(stan_file=stan_file)\n",
    "\n",
    "print(model.code())"
   ]
  },
  {
   "cell_type": "markdown",
   "metadata": {},
   "source": [
    "### Run Stan's variational inference algorithm, obtain fitted estimates\n",
    "\n",
    "The `CmdStanModel` method `variational` runs CmdStan's ADVI algorithm.\n",
    "Because this algorithm is unstable and may fail to converge, we run it with argument `require_converged` set to `False`.  We also specify a seed, to avoid instabilities as well as for reproducibility."
   ]
  },
  {
   "cell_type": "code",
   "execution_count": null,
   "metadata": {},
   "outputs": [],
   "source": [
    "vb_fit = model.variational(data=data_file, require_converged=False, seed=123)"
   ]
  },
  {
   "cell_type": "markdown",
   "metadata": {},
   "source": [
    "The ADVI algorithm provides estimates of all model parameters.\n",
    "\n",
    "The `variational` method returns a `CmdStanVB` object, with method `stan_variables`, which\n",
    "returns the approximate estimates of all model parameters as a Python dictionary."
   ]
  },
  {
   "cell_type": "code",
   "execution_count": null,
   "metadata": {},
   "outputs": [],
   "source": [
    "print(vb_fit.stan_variables())"
   ]
  },
  {
   "cell_type": "markdown",
   "metadata": {},
   "source": [
    "Posteriordb provides reference posteriors for all models.  For the blr model, conditioned on the dataset `sblri.json`, the reference posteriors are in file [sblri-blr.json](https://github.com/stan-dev/posteriordb/blob/master/posterior_database/reference_posteriors/summary_statistics/mean/mean/sblri-blr.json)\n",
    "\n",
    "The reference posteriors for all elements of `beta` and `sigma` are all very close to $1.0$.\n",
    "\n",
    "The experiments reported in the paper [Pathfinder: Parallel quasi-Newton variational inference](https://arxiv.org/abs/2108.03782) by Zhang et al. show that mean-field ADVI provides a better estimate of the posterior, as measured by the 1-Wasserstein distance to the reference posterior, than 75 iterations of the warmup Phase I algorithm used by the NUTS-HMC sampler, furthermore, ADVI is more computationally efficient, requiring fewer evaluations of the log density and gradient functions.  Therefore, using the estimates from ADVI to initialize the parameter values for the NUTS-HMC sampler will allow the sampler to do a better job of adapting the stepsize and metric during warmup, resulting in better performance and estimation.\n"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": null,
   "metadata": {},
   "outputs": [],
   "source": [
    "vb_vars = vb_fit.stan_variables()\n",
    "mcmc_vb_inits_fit = model.sample(\n",
    "    data=data_file, inits=vb_vars, iter_warmup=75, seed=12345\n",
    ")"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": null,
   "metadata": {},
   "outputs": [],
   "source": [
    "mcmc_vb_inits_fit.summary()"
   ]
  },
  {
   "cell_type": "markdown",
   "metadata": {},
   "source": [
    "The sampler estimates match the reference posterior."
   ]
  },
  {
   "cell_type": "code",
   "execution_count": null,
   "metadata": {},
   "outputs": [],
   "source": [
    "print(mcmc_vb_inits_fit.diagnose())"
   ]
  },
  {
   "cell_type": "markdown",
   "metadata": {},
   "source": [
    "Using the default random parameter initializations, we need to run more warmup iteratons. If we only run 75 warmup iterations with random inits, the result fails to estimate `sigma` correctly.  It is necessary to run the model with at least 150 warmup iterations to produce a good set of estimates."
   ]
  },
  {
   "cell_type": "code",
   "execution_count": null,
   "metadata": {},
   "outputs": [],
   "source": [
    "mcmc_random_inits_fit = model.sample(data=data_file, iter_warmup=75, seed=12345)"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": null,
   "metadata": {},
   "outputs": [],
   "source": [
    "mcmc_random_inits_fit.summary()"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": null,
   "metadata": {},
   "outputs": [],
   "source": [
    "print(mcmc_random_inits_fit.diagnose())"
   ]
  }
 ],
 "metadata": {
  "kernelspec": {
   "display_name": "Python 3",
   "language": "python",
   "name": "python3"
  },
  "language_info": {
   "codemirror_mode": {
    "name": "ipython",
    "version": 3
   },
   "file_extension": ".py",
   "mimetype": "text/x-python",
   "name": "python",
   "nbconvert_exporter": "python",
   "pygments_lexer": "ipython3",
   "version": "3.8.5"
  }
 },
 "nbformat": 4,
 "nbformat_minor": 4
}
